Introduction

This notebook demostrates the core functionality of pymatgen, including the core objects representing Elements, Species, Lattices, and Structures.

Written using:

  • pymatgen==2018.3.13

By convention, we import pymatgen as mg.


In [1]:
import pymatgen as mg

Basic Element, Specie and Composition objects

Pymatgen contains a set of core classes to represent an Element, Specie and Composition. These objects contains useful properties such as atomic mass, ionic radii, etc. These core classes are loaded by default with pymatgen. An Element can be created as follows:


In [2]:
si = mg.Element("Si")
print("Atomic mass of Si is {}".format(si.atomic_mass))
print("Si has a melting point of {}".format(si.melting_point))
print("Ionic radii for Si: {}".format(si.ionic_radii))


Atomic mass of Si is 28.0855 amu
Si has a melting point of 1687.0 K
Ionic radii for Si: {4: 0.54}

You can see that units are printed for atomic masses and ionic radii. Pymatgen comes with a complete system of managing units in pymatgen.core.unit. A Unit is a subclass of float that attaches units and handles conversions. For example,


In [3]:
print("Atomic mass of Si in kg: {}".format(si.atomic_mass.to("kg")))


Atomic mass of Si in kg: 4.6637069207919995e-26 kg

Please refer to the Units example for more information on units. Species are like Elements, except they have an explicit oxidation state. They can be used wherever Element is used for the most part.


In [4]:
fe2 = mg.Specie("Fe", 2)
print(fe2.atomic_mass)
print(fe2.ionic_radius)


55.845 amu
0.92 ang

A Composition is essentially an immutable mapping of Elements/Species with amounts, and useful properties like molecular weight, get_atomic_fraction, etc. Note that you can conveniently either use an Element/Specie object or a string as keys (this is a feature).


In [5]:
comp = mg.Composition("Fe2O3")
print("Weight of Fe2O3 is {}".format(comp.weight))
print("Amount of Fe in Fe2O3 is {}".format(comp["Fe"]))
print("Atomic fraction of Fe is {}".format(comp.get_atomic_fraction("Fe")))
print("Weight fraction of Fe is {}".format(comp.get_wt_fraction("Fe")))


Weight of Fe2O3 is 159.6882 amu
Amount of Fe in Fe2O3 is 2.0
Atomic fraction of Fe is 0.4
Weight fraction of Fe is 0.699425505453753 

Lattice & Structure objects

A Lattice represents a Bravais lattice. Convenience static functions are provided for the creation of common lattice types from a minimum number of arguments.


In [6]:
# Creates cubic Lattice with lattice parameter 4.2
lattice = mg.Lattice.cubic(4.2)
print(lattice.lengths_and_angles)


((4.2000000000000002, 4.2000000000000002, 4.2000000000000002), (90.0, 90.0, 90.0))

A Structure object represents a crystal structure (lattice + basis). A Structure is essentially a list of PeriodicSites with the same Lattice. Let us now create a CsCl structure.


In [7]:
structure = mg.Structure(lattice, ["Cs", "Cl"], [[0, 0, 0], [0.5, 0.5, 0.5]])
print("Unit cell vol = {}".format(structure.volume))
print("First site of the structure is {}".format(structure[0]))


Unit cell vol = 74.08800000000001
First site of the structure is [ 0.  0.  0.] Cs

The Structure object contains many useful manipulation functions. Since Structure is essentially a list, it contains a simple pythonic API for manipulation its sites. Some examples are given below. Please note that there is an immutable version of Structure known as IStructure, for the use case where you really need to enforce that the structure does not change. Conversion between these forms of Structure can be performed using from_sites().


In [8]:
structure.make_supercell([2, 2, 1]) #Make a 3 x 2 x 1 supercell of the structure
del structure[0] #Remove the first site
structure.append("Na", [0,0,0]) #Append a Na atom.
structure[-1] = "Li" #Change the last added atom to Li.
structure[0] = "Cs", [0.01, 0.5, 0] #Shift the first atom by 0.01 in fractional coordinates in the x-direction.
immutable_structure = mg.IStructure.from_sites(structure) #Create an immutable structure (cannot be modified).
print(immutable_structure)


Full Formula (Cs3 Li1 Cl4)
Reduced Formula: Cs3LiCl4
abc   :   8.400000   8.400000   4.200000
angles:  90.000000  90.000000  90.000000
Sites (8)
  #  SP       a     b    c
---  ----  ----  ----  ---
  0  Cs    0.01  0.5   0
  1  Cs    0.5   0     0
  2  Cs    0.5   0.5   0
  3  Cl    0.25  0.25  0.5
  4  Cl    0.25  0.75  0.5
  5  Cl    0.75  0.25  0.5
  6  Cl    0.75  0.75  0.5
  7  Li    0     0     0

Basic analyses

Pymatgen provides many analyses functions for Structures. Some common ones are given below.


In [9]:
#Determining the symmetry
from pymatgen.symmetry.analyzer import SpacegroupAnalyzer
finder = SpacegroupAnalyzer(structure)
print("The spacegroup is {}".format(finder.get_space_group_symbol()))


The spacegroup is Pmm2

We also have an extremely powerful structure matching tool.


In [10]:
from pymatgen.analysis.structure_matcher import StructureMatcher
#Let's create two structures which are the same topologically, but with different elements, and one lattice is larger.
s1 = mg.Structure(lattice, ["Cs", "Cl"], [[0, 0, 0], [0.5, 0.5, 0.5]])
s2 = mg.Structure(mg.Lattice.cubic(5), ["Rb", "F"], [[0, 0, 0], [0.5, 0.5, 0.5]])
m = StructureMatcher()
print(m.fit_anonymous(s1, s2)) #Returns a mapping which maps s1 and s2 onto each other. Strict element fitting is also available.


True

Input/output

Pymatgen also provides IO support for various file formats in the pymatgen.io package. A convenient set of read_structure and write_structure functions are also provided which auto-detects several well-known formats.


In [11]:
#Convenient IO to various formats. Format is intelligently determined from file name and extension.
structure.to(filename="POSCAR")
structure.to(filename="CsCl.cif")

#Or if you just supply fmt, you simply get a string.
print(structure.to(fmt="poscar"))
print(structure.to(fmt="cif"))


Cs3 Li1 Cl4
1.0
8.400000 0.000000 0.000000
0.000000 8.400000 0.000000
0.000000 0.000000 4.200000
Cs Cl Li
3 4 1
direct
0.010000 0.500000 0.000000 Cs
0.500000 0.000000 0.000000 Cs
0.500000 0.500000 0.000000 Cs
0.250000 0.250000 0.500000 Cl
0.250000 0.750000 0.500000 Cl
0.750000 0.250000 0.500000 Cl
0.750000 0.750000 0.500000 Cl
0.000000 0.000000 0.000000 Li

# generated using pymatgen
data_Cs3LiCl4
_symmetry_space_group_name_H-M   'P 1'
_cell_length_a   8.40000000
_cell_length_b   8.40000000
_cell_length_c   4.20000000
_cell_angle_alpha   90.00000000
_cell_angle_beta   90.00000000
_cell_angle_gamma   90.00000000
_symmetry_Int_Tables_number   1
_chemical_formula_structural   Cs3LiCl4
_chemical_formula_sum   'Cs3 Li1 Cl4'
_cell_volume   296.352
_cell_formula_units_Z   1
loop_
 _symmetry_equiv_pos_site_id
 _symmetry_equiv_pos_as_xyz
  1  'x, y, z'
loop_
 _atom_site_type_symbol
 _atom_site_label
 _atom_site_symmetry_multiplicity
 _atom_site_fract_x
 _atom_site_fract_y
 _atom_site_fract_z
 _atom_site_occupancy
  Cs  Cs1  1  0.010000  0.500000  0.000000  1
  Cs  Cs2  1  0.500000  0.000000  0.000000  1
  Cs  Cs3  1  0.500000  0.500000  0.000000  1
  Cl  Cl4  1  0.250000  0.250000  0.500000  1
  Cl  Cl5  1  0.250000  0.750000  0.500000  1
  Cl  Cl6  1  0.750000  0.250000  0.500000  1
  Cl  Cl7  1  0.750000  0.750000  0.500000  1
  Li  Li8  1  0.000000  0.000000  0.000000  1


In [12]:
#Reading a structure from a file.
structure = mg.Structure.from_file("POSCAR")

The vaspio_set module provides a means o obtain a complete set of VASP input files for performing calculations. Several useful presets based on the parameters used in the Materials Project are provided.


In [13]:
from pymatgen.io.vasp.sets import MPRelaxSet
v = MPRelaxSet(structure)
v.write_input("MyInputFiles") #Writes a complete set of input files for structure to the directory MyInputFiles

This concludes this pymatgen tutorial. Please explore the usage pages on pymatgen.org for more information.